The smuthi.fields package

fields

This subpackage contains functionality that has to do with the representation of electromagnetic fields in spherical or plane vector wave functions. The __init__ module contains some helper functions (e.g. with respect to SVWF indexing) and is the place to store default coordinate arrays for Sommerfeld integrals and field expansions.

smuthi.fields.angular_arrays(angular_resolution=0.008726646259971648)

Return azimuthal and polar angular arrays with a certain angular resolution

smuthi.fields.angular_frequency(vacuum_wavelength)

Angular frequency \(\omega = 2\pi c / \lambda\)

Parameters:vacuum_wavelength (float) – Vacuum wavelength in length unit
Returns:Angular frequency in the units of c=1 (time units=length units). This is at the same time the vacuum wavenumber.
smuthi.fields.blocksize

Number of coefficients in outgoing or regular spherical wave expansion for a single particle.

Parameters:
  • l_max (int) – Maximal multipole degree
  • m_max (int) – Maximal multipole order
Returns:

Number of indices for one particle, which is the maximal index plus 1.

smuthi.fields.branchpoint_correction(layer_refractive_indices, n_effective_array, neff_minimal_branchpoint_distance)

Check if an array of complex effective refractive index values (e.g. for Sommerfeld integration) contains possible branchpoint singularities and if so, replace them by nearby non-singular locations.

Parameters:
  • layer_refractive_indices (list or array) – Complex refractive indices of planarly layered medium
  • n_effective_array (1d numpy.array) – Complex effective refractive indexc values that are to be checked for branchpoint collision This array is changed during the function evaluation!
  • neff_minimal_branchpoint_distance (float) – Minimal distance that contour points shall have from branchpoint singularities
Returns:

corrected n_effective_array

smuthi.fields.create_k_parallel_array(vacuum_wavelength, neff_waypoints, neff_resolution)

Construct an array of complex in-plane wavenumbers (i.e., the radial component of the cylindrical coordinates of the wave-vector). This is used for the plane wave expansion of fields and for Sommerfeld integrals. Complex contours are used to improve numerical stability (see section 3.10.2.1 of [Egel 2018 dissertation]).

Parameters:
  • vacuum_wavelength (float) – Vacuum wavelength \(\lambda\) (length)
  • neff_waypoints (list or ndarray) – Corner points through which the contour runs This quantity is dimensionless (effective refractive index, will be multiplied by vacuum wavenumber)
  • neff_resolution (float) – Resolution of contour, again in terms of effective refractive index
Returns:

Array \(\kappa_i\) of in-plane wavenumbers (inverse length)

smuthi.fields.create_neff_array(neff_waypoints, neff_resolution)

Construct an array of complex effective refractive index values. The effective refractive index is a dimensionless quantity that will be multiplied by vacuum wavenumber to yield the in-plane component of a wave vector. This is used for the plane wave expansion of fields and for Sommerfeld integrals. Complex contours are used to improve numerical stability (see section 3.10.2.1 of [Egel 2018 dissertation]).

Parameters:
  • neff_waypoints (list or ndarray) – Corner points through which the contour runs
  • neff_resolution (float) – Resolution of contour (i.e., distance between adjacent elements)
Returns:

Array of complex effective refractive index values

smuthi.fields.default_Sommerfeld_k_parallel_array = None

Default n_effective array for the initial field (beams, dipoles) - needs to be set, e.g. at beginning of simulation

smuthi.fields.default_polar_angles = None

Default n_effective array for Sommerfeld integrals - needs to be set, e.g. at beginning of simulation

smuthi.fields.k_z(k_parallel=None, n_effective=None, k=None, omega=None, vacuum_wavelength=None, refractive_index=None)

z-component \(k_z=\sqrt{k^2-\kappa^2}\) of the wavevector. The branch cut is defined such that the imaginary part is not negative, compare section 2.3.1 of [Egel 2018 dissertation]. Not all of the arguments need to be specified.

Parameters:
  • k_parallel (numpy ndarray) – In-plane wavenumber \(\kappa\) (inverse length)
  • n_effective (numpy ndarray) – Effective refractive index \(n_\mathrm{eff}\)
  • k (float) – Wavenumber (inverse length)
  • omega (float) – Angular frequency \(\omega\) or vacuum wavenumber (inverse length, c=1)
  • vacuum_wavelength (float) – Vacuum wavelength \(\lambda\) (length)
  • refractive_index (complex) – Refractive index \(n_i\) of material
Returns:

z-component \(k_z\) of wavenumber with non-negative imaginary part (inverse length)

smuthi.fields.multi_to_single_index

Unique single index for the totality of indices characterizing a svwf expansion coefficient.

The mapping follows the scheme:

single index spherical wave expansion indices
\(n\) \(\tau\) \(l\) \(m\)
1 1 1 -1
2 1 1 0
3 1 1 1
4 1 2 -2
5 1 2 -1
6 1 2 0
1 l_max m_max
2 1 -1
Parameters:
  • tau (int) – Polarization index \(\tau\) (0=spherical TE, 1=spherical TM)
  • l (int) – Degree \(l\) (1, …, lmax)
  • m (int) – Order \(m\) (-min(l,mmax),…,min(l,mmax))
  • l_max (int) – Maximal multipole degree
  • m_max (int) – Maximal multipole order
Returns:

single index (int) subsuming \((\tau, l, m)\)

smuthi.fields.reasonable_Sommerfeld_kpar_contour(vacuum_wavelength, neff_waypoints=None, layer_refractive_indices=None, neff_imag=0.01, neff_max=None, neff_max_offset=1, neff_resolution=0.01, neff_minimal_branchpoint_distance=None)

Return a reasonable k_parallel array that is suitable as a Sommerfeld integral contour. Use this function if you don’t want to care for numerical details of your simulation.

Parameters:
  • vacuum_wavelength (float) – Vacuum wavelength \(\lambda\) (length)
  • neff_waypoints (list or ndarray) – Corner points through which the contour runs This quantity is dimensionless (effective refractive index, will be multiplied by vacuum wavenumber) If not provided, reasonable waypoints are estimated.
  • layer_refractive_indices (list) – Complex refractive indices of planarly layered medium Only needed when no neff_waypoints are provided
  • neff_imag (float) – Extent of the contour into the negative imaginary direction (in terms of effective refractive index, n_eff=kappa/omega). Only needed when no neff_waypoints are provided
  • neff_max (float) – Truncation value of contour (in terms of effective refractive index). Only needed when no neff_waypoints are provided
  • neff_max_offset (float) – Use the last estimated singularity location plus this value (in terms of effective refractive index). Default=1 Only needed when no neff_waypoints are provided and if no value for neff_max is specified.
  • neff_resolution (float) – Resolution of contour, again in terms of effective refractive index
  • neff_minimal_branchpoint_distance (float) – Minimal distance that contour points shall have from branchpoint singularities (in terms of effective refractive index). This is only relevant if not deflected into imaginary. Default: One fifth of neff_resolution
Returns:

Array \(\kappa_i\) of in-plane wavenumbers (inverse length)

smuthi.fields.reasonable_Sommerfeld_neff_contour(neff_waypoints=None, layer_refractive_indices=None, neff_imag=0.01, neff_max=None, neff_max_offset=1, neff_resolution=0.01, neff_minimal_branchpoint_distance=None)

Return a reasonable n_effective array that is suitable for the construction of a Sommerfeld k_parallel integral contour. Use this function if you don’t want to care for numerical details of your simulation.

Parameters:
  • neff_waypoints (list or ndarray) – Corner points through which the contour runs This quantity is dimensionless (effective refractive index, will be multiplied by vacuum wavenumber) If not provided, reasonable waypoints are estimated.
  • layer_refractive_indices (list) – Complex refractive indices of planarly layered medium Only needed when no neff_waypoints are provided
  • neff_imag (float) – Extent of the contour into the negative imaginary direction (in terms of effective refractive index, n_eff=kappa/omega). Only needed when no neff_waypoints are provided
  • neff_max (float) – Truncation value of contour (in terms of effective refractive index). Only needed when no neff_waypoints are provided
  • neff_max_offset (float) – Use the last estimated singularity location plus this value (in terms of effective refractive index). Default=1 Only needed when no neff_waypoints are provided and if no value for neff_max is specified.
  • neff_resolution (float) – Resolution of contour, again in terms of effective refractive index
  • neff_minimal_branchpoint_distance (float) – Minimal distance that contour points shall have from branchpoint singularities (in terms of effective refractive index). This is only relevant if not deflected into imaginary. Default: One fifth of neff_resolution
Returns:

Array of complex effective refractive index values

smuthi.fields.reasonable_neff_waypoints(layer_refractive_indices=None, neff_imag=0.01, neff_max=None, neff_max_offset=1)

Construct a reasonable list of waypoints for a k_parallel array of plane wave expansions. The waypoints mark a contour through the complex plane such that possible waveguide mode and branchpoint singularity locations are avoided (see section 3.10.2.1 of [Egel 2018 dissertation]).

Parameters:
  • layer_refractive_indices (list or array) – Complex refractive indices of the plane layer system
  • neff_imag (float) – Extent of the contour into the negative imaginary direction (in terms of effective refractive index, n_eff=kappa/omega).
  • neff_max (float) – Truncation value of contour (in terms of effective refractive index).
  • neff_max_offset (float) – If no value for neff_max is specified, use the last estimated singularity location plus this value (in terms of effective refractive index). Default=1
Returns:

List of complex waypoint values.

fields.expansions

Classes to manage the expansion of the electric field in plane wave and spherical wave basis sets.

class smuthi.fields.expansions.FieldExpansion

Base class for field expansions.

diverging(x, y, z)

Test if points are in domain where expansion could diverge. Virtual method to be overwritten in child classes.

Parameters:
  • x (numpy.ndarray) – x-coordinates of query points
  • y (numpy.ndarray) – y-coordinates of query points
  • z (numpy.ndarray) – z-coordinates of query points
Returns:

numpy.ndarray of bool datatype indicating if points are inside divergence domain.

electric_field(x, y, z)

Evaluate electric field. Virtual method to be overwritten in child classes.

Parameters:
  • x (numpy.ndarray) – x-coordinates of query points
  • y (numpy.ndarray) – y-coordinates of query points
  • z (numpy.ndarray) – z-coordinates of query points
Returns:

Tuple of (E_x, E_y, E_z) numpy.ndarray objects with the Cartesian coordinates of complex electric field.

magnetic_field(x, y, z, vacuum_wavelength)

Evaluate magnetic field. Virtual method to be overwritten in child classes.

Parameters:
  • x (numpy.ndarray) – x-coordinates of query points
  • y (numpy.ndarray) – y-coordinates of query points
  • z (numpy.ndarray) – z-coordinates of query points
  • vacuum_wavelength (float) – Vacuum wavelength in length units
Returns:

Tuple of (H_x, H_y, H_z) numpy.ndarray objects with the Cartesian coordinates of complex magnetic field.

valid(x, y, z)

Test if points are in definition range of the expansion. Abstract method to be overwritten in child classes.

Parameters:
  • x (numpy.ndarray) – x-coordinates of query points
  • y (numpy.ndarray) – y-coordinates of query points
  • z (numpy.ndarray) – z-coordinates of query points
Returns:

numpy.ndarray of bool datatype indicating if points are inside definition domain.

class smuthi.fields.expansions.PiecewiseFieldExpansion

Manage a field that is expanded in different ways for different domains, i.e., an expansion of the kind

\[\mathbf{E}(\mathbf{r}) = \sum_{i} \mathbf{E}_i(\mathbf{r}),\]

where

\[\begin{split}\mathbf{E}_i(\mathbf{r}) = \begin{cases} \tilde{\mathbf{E}}_i(\mathbf{r}) & \text{ if }\mathbf{r}\in D_i \\ 0 & \text{ else} \end{cases}\end{split}\]

and \(\tilde{\mathbf{E_i}}(\mathbf{r})\) is either a plane wave expansion or a spherical wave expansion, and \(D_i\) is its domain of validity.

compatible(other)

Returns always true, because any field expansion can be added to a piecewise field expansion.

diverging(x, y, z)

Test if points are in domain where expansion could diverge.

Parameters:
  • x (numpy.ndarray) – x-coordinates of query points
  • y (numpy.ndarray) – y-coordinates of query points
  • z (numpy.ndarray) – z-coordinates of query points
Returns:

numpy.ndarray of bool datatype indicating if points are inside divergence domain.

electric_field(x, y, z)

Evaluate electric field.

Parameters:
  • x (numpy.ndarray) – x-coordinates of query points
  • y (numpy.ndarray) – y-coordinates of query points
  • z (numpy.ndarray) – z-coordinates of query points
Returns:

Tuple of (E_x, E_y, E_z) numpy.ndarray objects with the Cartesian coordinates of complex electric field.

magnetic_field(x, y, z, vacuum_wavelength)

Evaluate magnetic field. Virtual method to be overwritten in child classes.

Parameters:
  • x (numpy.ndarray) – x-coordinates of query points
  • y (numpy.ndarray) – y-coordinates of query points
  • z (numpy.ndarray) – z-coordinates of query points
  • vacuum_wavelength (float) – Vacuum wavelength in length units
Returns:

Tuple of (H_x, H_y, H_z) numpy.ndarray objects with the Cartesian coordinates of complex magnetic field.

valid(x, y, z)

Test if points are in definition range of the expansion.

Parameters:
  • x (numpy.ndarray) – x-coordinates of query points
  • y (numpy.ndarray) – y-coordinates of query points
  • z (numpy.ndarray) – z-coordinates of query points
Returns:

numpy.ndarray of bool datatype indicating if points are inside definition domain.

class smuthi.fields.expansions.PlaneWaveExpansion(k, k_parallel, azimuthal_angles, kind=None, reference_point=None, lower_z=-inf, upper_z=inf)

A class to manage plane wave expansions of the form

\[\mathbf{E}(\mathbf{r}) = \sum_{j=1}^2 \iint \mathrm{d}^2\mathbf{k}_\parallel \, g_{j}(\kappa, \alpha) \mathbf{\Phi}^\pm_j(\kappa, \alpha; \mathbf{r} - \mathbf{r}_i)\]

for \(\mathbf{r}\) located in a layer defined by \(z\in [z_{min}, z_{max}]\) and \(\mathrm{d}^2\mathbf{k}_\parallel = \kappa\,\mathrm{d}\alpha\,\mathrm{d}\kappa\).

The double integral runs over \(\alpha\in[0, 2\pi]\) and \(\kappa\in[0,\kappa_\mathrm{max}]\). Further, \(\mathbf{\Phi}^\pm_j\) are the PVWFs, see plane_vector_wave_function().

Internally, the expansion coefficients \(g_{ij}^\pm(\kappa, \alpha)\) are stored as a 3-dimensional numpy ndarray.

If the attributes k_parallel and azimuthal_angles have only a single entry, a discrete distribution is assumed:

\[g_{j}^\pm(\kappa, \alpha) \sim \delta^2(\mathbf{k}_\parallel - \mathbf{k}_{\parallel, 0})\]
Parameters:
  • k (float) – wavenumber in layer where expansion is valid
  • k_parallel (numpy ndarray) – array of in-plane wavenumbers (can be float or complex)
  • azimuthal_angles (numpy ndarray) – \(\alpha\), from 0 to \(2\pi\)
  • kind (str) – ‘upgoing’ for \(g^+\) and ‘downgoing’ for \(g^-\) type expansions
  • reference_point (list or tuple) – [x, y, z]-coordinates of point relative to which the plane waves are defined.
  • lower_z (float) – the expansion is valid on and above that z-coordinate
  • upper_z (float) – the expansion is valid below that z-coordinate
coefficients

coefficients[j, k, l] contains

Type:numpy ndarray
:math:`g^pm_{j}
Type:kappa_{k}, alpha_{l}
class OptimizationMethodsForLinux

An enumeration.

evaluate_r_times_eikr

Attention! Sometimes this function can decrease speed on 1 core mode. Here foo_x, foo_y, foo_z are supposed to be 2dim arrays with [None, :, :]. This function can replace snippet

exp_j = np.exp(1j * exp_feed) foo_x_eikr = foo_x * exp_j foo_y_eikr = foo_y * exp_j foo_z_eikr = foo_z * exp_j

by

foo_x_eikr, foo_y_eikr, foo_z_eikr = numba_multiple_on_exp(foo_x, foo_y, foo_z, kr).

numba_3tensordots_1dim_times_2dim

This function can replace snippet

‘foo = np.tensordot(x_float_1dim, x_complex_2dim, axes=0)

foo += np.tensordot(y_float_1dim, y_complex_2dim, axes=0)

foo += np.tensordot(z_float_1dim, z_complex_2dim, axes=0)’

by

‘foo = get_3_tensordots(x_float_1dim, y_float_1dim, z_float_1dim,
x_complex_2dim, y_complex_2dim, z_complex_2dim)’
numba_trapz_3dim_array

This function can replace snippet

‘foo = np.trapz(y, x)’

by

‘foo = numba_trapz_3dim_array(y, x)’

class OptimizationMethodsFor_Not_Linux

An enumeration.

evaluate_r_times_eikr
numba_3tensordots_1dim_times_2dim
numba_trapz_3dim_array
class RawSliceOfField(axis, chunks, values)
azimuthal_angle_grid()

Meshgrid of azimuthal_angles with respect to n_effective

compatible(other)

Check if two plane wave expansions are compatible in the sense that they can be added coefficient-wise

Parameters:other (FieldExpansion) – expansion object to add to this object
Returns:bool (true if compatible, false else)
diverging(x, y, z)

Test if points are in domain where expansion could diverge.

Parameters:
  • x (numpy.ndarray) – x-coordinates of query points
  • y (numpy.ndarray) – y-coordinates of query points
  • z (numpy.ndarray) – z-coordinates of query points
Returns:

numpy.ndarray of bool datatype indicating if points are inside divergence domain.

electric_field(x, y, z, max_chunksize=50, cpu_precision='single precision')

Evaluate electric field.

Parameters:
  • x (numpy.ndarray) – x-coordinates of query points
  • y (numpy.ndarray) – y-coordinates of query points
  • z (numpy.ndarray) – z-coordinates of query points
  • max_chunksize (int) – max number of field points that are simultaneously evaluated when running in CPU mode. In Windows/MacOS max_chunksize = chunksize, in Linux it can be decreased considering available CPU cores.
  • cpu_precision (string) – set ‘double precision’ to use float64 and complex128 types instead of float32 and complex64
Returns:

Tuple of (E_x, E_y, E_z) numpy.ndarray objects with the Cartesian coordinates of complex electric field.

k_parallel_grid()

Meshgrid of n_effective with respect to azimuthal_angles

k_z()
k_z_grid()
magnetic_field(x, y, z, vacuum_wavelength, max_chunksize=50, cpu_precision='single precision')

Evaluate magnetic field.

Parameters:
  • x (numpy.ndarray) – x-coordinates of query points
  • y (numpy.ndarray) – y-coordinates of query points
  • z (numpy.ndarray) – z-coordinates of query points
  • vacuum_wavelength (float) – Vacuum wavelength in length units
  • chunksize (int) – number of field points that are simultaneously evaluated when running in CPU mode
Returns:

Tuple of (H_x, H_y, H_z) numpy.ndarray objects with the Cartesian coordinates of complex magnetic field.

set_reference_point(new_reference_point)

Set a new reference point. This implies also a phase factor on the coefficients.

Parameters:new_reference_point (list or tuple) – [x, y, z]-coordinates of point relative to which the plane waves are defined.
valid(x, y, z)

Test if points are in definition range of the expansion.

Parameters:
  • x (numpy.ndarray) – x-coordinates of query points
  • y (numpy.ndarray) – y-coordinates of query points
  • z (numpy.ndarray) – z-coordinates of query points
Returns:

numpy.ndarray of bool datatype indicating if points are inside definition domain.

class smuthi.fields.expansions.SphericalWaveExpansion(k, l_max, m_max=None, kind=None, reference_point=None, lower_z=-inf, upper_z=inf, inner_r=0, outer_r=inf)

A class to manage spherical wave expansions of the form

\[\mathbf{E}(\mathbf{r}) = \sum_{\tau=1}^2 \sum_{l=1}^\infty \sum_{m=-l}^l a_{\tau l m} \mathbf{\Psi}^{(\nu)}_{\tau l m}(\mathbf{r} - \mathbf{r}_i)\]

for \(\mathbf{r}\) located in a layer defined by \(z\in [z_{min}, z_{max}]\) and where \(\mathbf{\Psi}^{(\nu)}_{\tau l m}\) are the SVWFs, see smuthi.vector_wave_functions.spherical_vector_wave_function().

Internally, the expansion coefficients \(a_{\tau l m}\) are stored as a 1-dimensional array running over a multi index \(n\) subsumming over the SVWF indices \((\tau,l,m)\). The mapping from the SVWF indices to the multi index is organized by the function multi_to_single_index().

Parameters:
  • k (float) – wavenumber in layer where expansion is valid
  • l_max (int) – maximal multipole degree \(l_\mathrm{max}\geq 1\) where to truncate the expansion.
  • m_max (int) – maximal multipole order \(0 \leq m_\mathrm{max} \leq l_\mathrm{max}\) where to truncate the expansion.
  • kind (str) – ‘regular’ for \(\nu=1\) or ‘outgoing’ for \(\nu=3\)
  • reference_point (list or tuple) – [x, y, z]-coordinates of point relative to which the spherical waves are considered (e.g., particle center).
  • lower_z (float) – the expansion is valid on and above that z-coordinate
  • upper_z (float) – the expansion is valid below that z-coordinate
  • inner_r (float) – radius inside which the expansion diverges (e.g. circumscribing sphere of particle)
  • outer_r (float) – radius outside which the expansion diverges
coefficients

expansion coefficients \(a_{\tau l m}\) ordered by multi index \(n\)

Type:numpy ndarray
coefficients_tlm(tau, l, m)

SWE coefficient for given (tau, l, m)

Parameters:
  • tau (int) – SVWF polarization (0 for spherical TE, 1 for spherical TM)
  • l (int) – SVWF degree
  • m (int) – SVWF order
Returns:

SWE coefficient

compatible(other)

Check if two spherical wave expansions are compatible in the sense that they can be added coefficient-wise

Parameters:other (FieldExpansion) – expansion object to add to this object
Returns:bool (true if compatible, false else)
diverging(x, y, z)

Test if points are in domain where expansion could diverge.

Parameters:
  • x (numpy.ndarray) – x-coordinates of query points
  • y (numpy.ndarray) – y-coordinates of query points
  • z (numpy.ndarray) – z-coordinates of query points
Returns:

numpy.ndarray of bool datatype indicating if points are inside divergence domain.

electric_field(x, y, z)

Evaluate electric field.

Parameters:
  • x (numpy.ndarray) – x-coordinates of query points
  • y (numpy.ndarray) – y-coordinates of query points
  • z (numpy.ndarray) – z-coordinates of query points
Returns:

Tuple of (E_x, E_y, E_z) numpy.ndarray objects with the Cartesian coordinates of complex electric field.

magnetic_field(x, y, z, vacuum_wavelength)

Evaluate magnetic field.

Parameters:
  • x (numpy.ndarray) – x-coordinates of query points
  • y (numpy.ndarray) – y-coordinates of query points
  • z (numpy.ndarray) – z-coordinates of query points
  • vacuum_wavelength (float) – Vacuum wavelength in length units
Returns:

Tuple of (H_x, H_y, H_z) numpy.ndarray objects with the Cartesian coordinates of complex electric field.

valid(x, y, z)

Test if points are in definition range of the expansion.

Parameters:
  • x (numpy.ndarray) – x-coordinates of query points
  • y (numpy.ndarray) – y-coordinates of query points
  • z (numpy.ndarray) – z-coordinates of query points
Returns:

numpy.ndarray of bool datatype indicating if points are inside definition domain.

fields.expansions_cuda

This module contains CUDA source code for the evaluation of the electric field from a VWF expansion.

fields.transformatinos

Functions for the transformation of plane and spherical vector wave functions as well as of plane and spherical wave fex.

smuthi.fields.transformations.block_rotation_matrix_D_svwf(l_max, m_max, alpha, beta, gamma, wdsympy=False)

Rotation matrix for the rotation of SVWFs between the labratory coordinate system (L) and a rotated coordinate system (R)

Parameters:
  • l_max (int) – Maximal multipole degree
  • m_max (int) – Maximal multipole order
  • alpha (float) – First Euler angle, rotation around z-axis, in rad
  • beta (float) – Second Euler angle, rotation around y’-axis in rad
  • gamma (float) – Third Euler angle, rotation around z’’-axis in rad
  • wdsympy (bool) – If True, Wigner-d-functions come from the sympy toolbox
Returns:

rotation matrix of dimension [blocksize, blocksize]

smuthi.fields.transformations.pwe_to_swe_conversion(pwe, l_max, m_max, reference_point)

Convert plane wave expansion object to a spherical wave expansion object.

Parameters:
  • pwe (PlaneWaveExpansion) – Plane wave expansion to be converted
  • l_max (int) – Maximal multipole degree of spherical wave expansion
  • m_max (int) – Maximal multipole order of spherical wave expansion
  • reference_point (list) – Coordinates of reference point in the format [x, y, z]
Returns:

SphericalWaveExpansion object.

smuthi.fields.transformations.swe_to_pwe_conversion(swe, k_parallel, azimuthal_angles, layer_system=None, layer_number=None, layer_system_mediated=False, only_l=None, only_m=None, only_pol=None, only_tau=None)

Convert SphericalWaveExpansion object to a PlaneWaveExpansion object.

Parameters:
  • swe (SphericalWaveExpansion) – Spherical wave expansion to be converted
  • k_parallel (numpy array or str) – In-plane wavenumbers for the pwe object.
  • azimuthal_angles (numpy array or str) – Azimuthal angles for the pwe object
  • layer_system (smuthi.layers.LayerSystem) – Stratified medium in which the origin of the SWE is located
  • layer_number (int) – Layer number in which the PWE should be valid.
  • layer_system_mediated (bool) – If True, the PWE refers to the layer system response of the SWE, otherwise it is the direct transform.
  • only_pol (int) – if set to 0 or 1, only this plane wave polarization (0=TE, 1=TM) is considered
  • only_tau (int) – if set to 0 or 1, only this spherical vector wave polarization (0 — magnetic, 1 — electric) is considered
  • only_l (int) – if set to positive number, only this multipole degree is considered
  • only_m (int) – if set to non-negative number, only this multipole order is considered
Returns:

Tuple of two PlaneWaveExpansion objects, first upgoing, second downgoing.

smuthi.fields.transformations.transformation_coefficients_vwf(tau, l, m, pol, kp=None, kz=None, pilm_list=None, taulm_list=None, dagger=False)

Transformation coefficients B to expand SVWF in PVWF and vice versa:

\[B_{\tau l m,j}(x) = -\frac{1}{\mathrm{i}^{l+1}} \frac{1}{\sqrt{2l(l+1)}} (\mathrm{i} \delta_{j1} + \delta_{j2}) (\delta_{\tau j} \tau_l^{|m|}(x) + (1-\delta_{\tau j} m \pi_l^{|m|}(x))\]

For the definition of the \(\tau_l^m\) and \(\pi_l^m\) functions, see A. Doicu, T. Wriedt, and Y. A. Eremin: “Light Scattering by Systems of Particles”, Springer-Verlag, 2006

Compare also section 2.3.3 of [Egel 2018 diss].

Parameters:
  • tau (int) – SVWF polarization, 0 for spherical TE, 1 for spherical TM
  • l (int) – l=1,… SVWF multipole degree
  • m (int) – m=-l,…,l SVWF multipole order
  • pol (int) – PVWF polarization, 0 for TE, 1 for TM
  • kp (numpy array) – PVWF in-plane wavenumbers
  • kz (numpy array) – complex numpy-array: PVWF out-of-plane wavenumbers
  • pilm_list (list) – 2D list numpy-arrays: alternatively to kp and kz, pilm and taulm as generated with legendre_normalized can directly be handed
  • taulm_list (list) – 2D list numpy-arrays: alternatively to kp and kz, pilm and taulm as generated with legendre_normalized can directly be handed
  • dagger (bool) – switch on when expanding PVWF in SVWF and off when expanding SVWF in PVWF
Returns:

Transformation coefficient as array (size like kp).

smuthi.fields.transformations.translation_block(vacuum_wavelength, receiving_particle, emitting_particle, layer_system, kind)
Direct particle translation matrix \(W\) for two particles that do not have intersecting circumscribing spheres.

This routine is explicit.

To reduce computation time, this routine relies on two internal accelerations. First, in most cases the number of unique maximum multipole indicies, \((\tau, l_{max}, m_{max})\), is much less than the number of unique particles. Therefore, all calculations that depend only on multipole indicies are stored in an intermediate hash table. Second, Cython acceleration is used by default to leverage fast looping. If the Cython files are not supported, this routine will fall back on equivalent Python looping.

Cython acceleration can be between 10-1,000x faster compared to the Python equivalent. Speed variability depends on the number of unique multipoles indicies, the size of the largest multipole order, and if particles share the same z coordinate.

Parameters:
Returns:

Direct coupling matrix block as numpy array.

smuthi.fields.transformations.translation_coefficients_svwf(tau1, l1, m1, tau2, l2, m2, k, d, sph_hankel=None, legendre=None, exp_immphi=None)

Coefficients of the translation operator for the expansion of an outgoing spherical wave in terms of regular spherical waves with respect to a different origin:

\[\mathbf{\Psi}_{\tau l m}^{(3)}(\mathbf{r} + \mathbf{d} = \sum_{\tau'} \sum_{l'} \sum_{m'} A_{\tau l m, \tau' l' m'} (\mathbf{d}) \mathbf{\Psi}_{\tau' l' m'}^{(1)}(\mathbf{r})\]

for \(|\mathbf{r}|<|\mathbf{d}|\).

See also section 2.3.3 and appendix B of [Egel 2018 diss].

Parameters:
  • tau1 (int) – tau1=0,1: Original wave’s spherical polarization
  • l1 (int) – l=1,…: Original wave’s SVWF multipole degree
  • m1 (int) – m=-l,…,l: Original wave’s SVWF multipole order
  • tau2 (int) – tau2=0,1: Partial wave’s spherical polarization
  • l2 (int) – l=1,…: Partial wave’s SVWF multipole degree
  • m2 (int) – m=-l,…,l: Partial wave’s SVWF multipole order
  • k (float or complex) – wavenumber (inverse length unit)
  • d (list) – translation vectors in format [dx, dy, dz] (length unit) dx, dy, dz can be scalars or ndarrays
  • sph_hankel (list) – Optional. sph_hankel[i] contains the spherical hankel funciton of degree i, evaluated at k*d where d is the norm of the distance vector(s)
  • legendre (list) – Optional. legendre[l][m] contains the legendre function of order l and degree m, evaluated at cos(theta) where theta is the polar angle(s) of the distance vector(s)
Returns:

translation coefficient A (complex)

smuthi.fields.transformations.translation_coefficients_svwf_out_to_out(tau1, l1, m1, tau2, l2, m2, k, d, sph_bessel=None, legendre=None, exp_immphi=None)

Coefficients of the translation operator for the expansion of an outgoing spherical wave in terms of outgoing spherical waves with respect to a different origin:

\[\mathbf{\Psi}_{\tau l m}^{(3)}(\mathbf{r} + \mathbf{d} = \sum_{\tau'} \sum_{l'} \sum_{m'} A_{\tau l m, \tau' l' m'} (\mathbf{d}) \mathbf{\Psi}_{\tau' l' m'}^{(3)}(\mathbf{r})\]

for \(|\mathbf{r}|>|\mathbf{d}|\).

Parameters:
  • tau1 (int) – tau1=0,1: Original wave’s spherical polarization
  • l1 (int) – l=1,…: Original wave’s SVWF multipole degree
  • m1 (int) – m=-l,…,l: Original wave’s SVWF multipole order
  • tau2 (int) – tau2=0,1: Partial wave’s spherical polarization
  • l2 (int) – l=1,…: Partial wave’s SVWF multipole degree
  • m2 (int) – m=-l,…,l: Partial wave’s SVWF multipole order
  • k (float or complex) – wavenumber (inverse length unit)
  • d (list) – translation vectors in format [dx, dy, dz] (length unit) dx, dy, dz can be scalars or ndarrays
  • sph_bessel (list) – Optional. sph_bessel[i] contains the spherical Bessel funciton of degree i, evaluated at k*d where d is the norm of the distance vector(s)
  • legendre (list) – Optional. legendre[l][m] contains the legendre function of order l and degree m, evaluated at cos(theta) where theta is the polar angle(s) of the distance vector(s)
Returns:

translation coefficient A (complex)

fields.vector_wave_functions

This module contains the vector wave functions and their transformations.

smuthi.fields.vector_wave_functions.plane_vector_wave_function(x, y, z, kp, alpha, kz, pol)

Electric field components of plane wave (PVWF), see section 2.3.1 of [Egel 2018 diss].

\[\mathbf{\Phi}_j = \exp ( \mathrm{i} \mathbf{k} \cdot \mathbf{r} ) \hat{ \mathbf{e} }_j\]

with \(\hat{\mathbf{e}}_0\) denoting the unit vector in azimuthal direction (‘TE’ or ‘s’ polarization), and \(\hat{\mathbf{e}}_1\) denoting the unit vector in polar direction (‘TM’ or ‘p’ polarization).

The input arrays should have one of the following dimensions:

  • x,y,z: (N x 1) matrix
  • kp,alpha,kz: (1 x M) matrix
  • Ex, Ey, Ez: (M x N) matrix

or

  • x,y,z: (M x N) matrix
  • kp,alpha,kz: scalar
  • Ex, Ey, Ez: (M x N) matrix
Parameters:
  • x (numpy.ndarray) – x-coordinate of position where to test the field (length unit)
  • y (numpy.ndarray) – y-coordinate of position where to test the field
  • z (numpy.ndarray) – z-coordinate of position where to test the field
  • kp (numpy.ndarray) – parallel component of k-vector (inverse length unit)
  • alpha (numpy.ndarray) – azimthal angle of k-vector (rad)
  • kz (numpy.ndarray) – z-component of k-vector (inverse length unit)
  • pol (int) – Polarization (0=TE, 1=TM)
Returns:

  • x-coordinate of PVWF electric field (numpy.ndarray)
  • y-coordinate of PVWF electric field (numpy.ndarray)
  • z-coordinate of PVWF electric field (numpy.ndarray)

smuthi.fields.vector_wave_functions.spherical_vector_wave_function(x, y, z, k, nu, tau, l, m)

Electric field components of spherical vector wave function (SVWF). The conventions are chosen according to A. Doicu, T. Wriedt, and Y. A. Eremin: “Light Scattering by Systems of Particles”, Springer-Verlag, 2006 See also section 2.3.2 of [Egel 2018 diss].

Parameters:
  • x (numpy.ndarray) – x-coordinate of position where to test the field (length unit)
  • y (numpy.ndarray) – y-coordinate of position where to test the field
  • z (numpy.ndarray) – z-coordinate of position where to test the field
  • k (float or complex) – wavenumber (inverse length unit)
  • nu (int) – 1 for regular waves, 3 for outgoing waves
  • tau (int) – spherical polarization, 0 for spherical TE and 1 for spherical TM
  • l (int) – l=1,… multipole degree (polar quantum number)
  • m (int) – m=-l,…,l multipole order (azimuthal quantum number)
Returns:

  • x-coordinate of SVWF electric field (numpy.ndarray)
  • y-coordinate of SVWF electric field (numpy.ndarray)
  • z-coordinate of SVWF electric field (numpy.ndarray)